(triv_lib/mrchcube.c:68)
Prototype:
MCPolygonStruct *MCThresholdCube(MCCubeCornerScalarStruct *CCS,
RealType Threshold)
Description:
Given 8 cube corner values (scalars), compute the polygon(s) in this cube
along the isosurface at Threshold. if CCS has gradient information, it is
used to approximate normals at the vertices.
7 K 6
***********************
* + * *
L * + * * Vertices 0 - 7
* + I * J * Edges A - L
4 *********************** 5 *
* + * *
* + * * G
* + H * *
* + * *
* + * F *
E * + C * *
* ++++++++++++++*+++++++* 2
* D + 3 * *
* + * * B
* + * *
***********************
0 A 1
Parameters:
CCS: | The cube's dimensions/information.
|
---|
Threshold: | Iso surface level.
|
---|
Returned Value:
MCPolygonStruct *: List of polygons (not necessarily triangles), or
possibly NULL.
|
---|
Keywords:
marching cubes
(triv_lib/triv_gen.c:93)
Prototype:
TrivTVStruct *TrivBspTVNew(int ULength,
int VLength,
int WLength,
int UOrder,
int VOrder,
int WOrder,
CagdPointType PType)
Description:
Allocates the memory required for a new Bspline trivariate.
Parameters:
ULength: | Number of control points in the U direction.
|
---|
VLength: | Number of control points in the V direction.
|
---|
WLength: | Number of control points in the W direction.
|
---|
UOrder: | Order of trivariate in the U direction.
|
---|
VOrder: | Order of trivariate in the V direction.
|
---|
WOrder: | Order of trivariate in the W direction.
|
---|
PType: | Type of control points (E2, P3, etc.).
|
---|
Returned Value:
TrivTVStruct *: An uninitialized freeform trivariate Bspline.
|
---|
Keywords:
trivariates
allocation
(triv_lib/trivrais.c:175)
Prototype:
TrivTVStruct *TrivBspTVDegreeRaise(TrivTVStruct *TV, TrivTVDirType Dir)
Description:
Returns a new Bspline surface, identical to the original but with one
degree higher, in the requested direction Dir.
Parameters:
TV: | To raise it degree by one.
|
---|
Dir: | Direction of degree raising. Either U, V or W.
|
---|
Returned Value:
TrivTVStruct *: A trivariate with one degree higher in direction Dir,
representing the same geometry as TV.
|
---|
Keywords:
degree raising
(triv_lib/triv_der.c:154)
Prototype:
TrivTVStruct *TrivBspTVDerive(TrivTVStruct *TV, TrivTVDirType Dir)
Description:
Given a Bspline trivariate, computes its partial derivative trivariate in
direction Dir.
Let old control polygon be P(i), i = 0 to k-1, and Q(i) be new one then:
Q(i) = (k - 1) * (P(i+1) - P(i)) / (Kv(i + k) - Kv(i + 1)), i = 0 to k-2.
Parameters:
TV: | Trivariate to differentiate.
|
---|
Dir: | Direction of differentiation. Either U or V or W.
|
---|
Returned Value:
TrivTVStruct *: Differentiated trivariate in direction Dir. A Bspline
trivariate.
|
---|
Keywords:
trivariates
(triv_lib/triv_ref.c:86)
Prototype:
TrivTVStruct *TrivBspTVKnotInsertNDiff(TrivTVStruct *TV,
TrivTVDirType Dir,
int Replace,
CagdRType *t,
int n)
Description:
Given a Bspline trivariate, inserts n knots with different values as
defined by t.
If, however, Replace is TRUE, the knot are simply replacing the current
knot vector in the prescribed direction.
Parameters:
TV: | Trivariate to refine according to t in direction Dir.
|
---|
Replace: | If TRUE t is a knot vector exaclt in the length of the knot
vector in direction Dir in TV and t simply replaces than knot
vector. If FALSE, the knot vector in direction Dir in TV is
refined by adding all the knots in t.
|
---|
t: | Knot vector to refine/replace the knot vector of TV in
direction Dir.
|
---|
n: | Length of vector t.
|
---|
Returned Value:
TrivTVStruct *: The refined trivariate. A Bspline trivariate.
|
---|
Keywords:
trivariates
(triv_lib/triv_gen.c:137)
Prototype:
TrivTVStruct *TrivBzrTVNew(int ULength,
int VLength,
int WLength,
CagdPointType PType)
Description:
Allocates the memory required for a new Bezier trivariate.
Parameters:
ULength: | Number of control points in the U direction.
|
---|
VLength: | Number of control points in the V direction.
|
---|
WLength: | Number of control points in the W direction.
|
---|
PType: | Type of control points (E2, P3, etc.).
|
---|
Returned Value:
TrivTVStruct *: An uninitialized freeform trivariate Bezier.
|
---|
Keywords:
trivariates
allocation
(triv_lib/trivrais.c:69)
Prototype:
TrivTVStruct *TrivBzrTVDegreeRaise(TrivTVStruct *TV, TrivTVDirType Dir)
Description:
Returns a new Bezier trivariate, identical to the original but with one
degree higher, in the requested direction Dir.
Let old control polygon be P(i), i = 0 to k-1, and Q(i) be new one then:
i k-i
Q(0) = P(0), Q(i) = --- P(i-1) + (---) P(i), Q(k) = P(k-1).
k k
This is applied to all rows/cols of the trivariate.
Parameters:
TV: | To raise it degree by one.
|
---|
Dir: | Direction of degree raising. Either U, V or W.
|
---|
Returned Value:
TrivTVStruct *: A surface with one degree higher in direction Dir,
representing the same geometry as TV.
|
---|
Keywords:
degree raising
(triv_lib/triv_der.c:66)
Prototype:
TrivTVStruct *TrivBzrTVDerive(TrivTVStruct *TV, TrivTVDirType Dir)
Description:
Given a Bezier trivariate, computes its partial derivative trivariate in
direction Dir.
Let old control polygon be P(i), i = 0 to k-1, and Q(i) be new one then:
Q(i) = (k - 1) * (P(i+1) - P(i)), i = 0 to k-2.
Parameters:
TV: | Trivariate to differentiate.
|
---|
Dir: | Direction of differentiation. Either U or V or W.
|
---|
Returned Value:
TrivTVStruct *: Differentiated trivariate in direction Dir. A Bezier
trivariate.
|
---|
Keywords:
trivariates
(triv_lib/triv_gen.c:507)
Prototype:
TrivTVStruct *TrivCnvrtBezier2BsplineTV(TrivTVStruct *TV)
Description:
Converts a Bezier trivariate into a Bspline trivariate by adding two open
end uniform knot vectors to it.
Parameters:
TV: | A Bezier trivariate to convert to a Bspline TV.
|
---|
Returned Value:
TrivTVStruct *: A Bspline trivariate representing the same geometry as
the given Bezier TV.
|
---|
Keywords:
conversion
trivariate
(triv_lib/trivcoer.c:25)
Prototype:
TrivTVStruct *TrivCoerceTVTo(TrivTVStruct *TV, CagdPointType PType)
Description:
Coerces a trivariate to point type PType.
Parameters:
TV: | To coerce to a new point type PType.
|
---|
PType: | New point type for TV.
|
---|
Returned Value:
TrivTVStruct *: A new trivariate with PType as its point type.
|
---|
Keywords:
coercion
(triv_lib/triv_dbg.c:25)
Prototype:
void TrivDbg(void *Obj)
Description:
Prints trivariates stderr. Should be linked to programs for debugging
purposes, so trivariates may be inspected from a debugger.
Parameters:
Obj: | A trivariate - to be printed to stderr.
|
---|
Returned Value:
Keywords:
debugging
(triv_lib/triv_err.c:59)
Prototype:
char *TrivDescribeError(TrivFatalErrorType ErrorNum)
Description:
Returns a string describing a the given error. Errors can be raised by
any member of this triv library as well as other users. Raised error will
cause an invokation of TrivFatalError function which decides how to handle
this error. TrivFatalError can for example, invoke this routine with the
error type, print the appropriate message and quit the program.
Parameters:
ErrorNum: | Type of the error that was raised.
|
---|
Returned Value:
char *: A string describing the error type.
|
---|
Keywords:
error handling
(triv_lib/triv_ftl.c:27)
Prototype:
void TrivFatalError(TrivFatalErrorType ErrID)
Description:
Trap Triv_lib errors right here. Provides a default error handler for the
triv library. Gets an error description using TrivDescribeError, prints it
and exit the program using exit.
Parameters:
ErrID: | Error type that was raised.
|
---|
Returned Value:
Keywords:
error handling
(triv_lib/trinterp.c:26)
Prototype:
TrivTVStruct *TrivInterpTrivar(TrivTVStruct *TV)
Description:
Interpolates control points of given trivariate, preserving the order
and continuity of the original trivariate.
Parameters:
TV: | Trivariate to interpolate its control mesh.
|
---|
Returned Value:
TrivTVStruct *: The interpolating trivariate.
|
---|
Keywords:
interpolation
(triv_lib/trivcmpt.c:50)
Prototype:
CagdBType TrivMakeTVsCompatible(TrivTVStruct **TV1,
TrivTVStruct **TV2,
CagdBType SameUOrder,
CagdBType SameVOrder,
CagdBType SameWOrder,
CagdBType SameUKV,
CagdBType SameVKV,
CagdBType SameWKV)
Description:
Given two trivariates, makes them compatible by:
1. Coercing their point type to be the same.
2. Making them have the same curve type.
3. Raising the degree of the lower one to be the same as the higher.
4. Refining them to a common knot vector (If Bspline and SameOrder).
Note 3 is performed if SameOrder TRUE, 4 if SameKV TRUE.
Both trivariates are modified IN PLACE.
Parameters:
TV1, TV2: | Two surfaces to be made compatible, in place.
|
---|
SameUOrder: | If TRUE, this routine make sure they share the same U
order.
|
---|
SameVOrder: | If TRUE, this routine make sure they share the same
order.
|
---|
SameWOrder: | If TRUE, this routine make sure they share the same W
order.
|
---|
SameUKV: | If TRUE, this routine make sure they share the same U
knot vector and hence continuity. *
|
---|
SameVKV: | If TRUE, this routine make sure they share the same
knot vector and hence continuity.
|
---|
SameWKV: | If TRUE, this routine make sure they share the same W
knot vector and hence continuity.
|
---|
Returned Value:
CagdBType: TRUE if successful, FALSE otherwise.
|
---|
Keywords:
compatibility
(triv_lib/triv_aux.c:82)
Prototype:
CagdBType TrivParamInDomain(TrivTVStruct *TV, CagdRType t, TrivTVDirType Dir)
Description:
Given a tri-variate and a domain - validate it.
Parameters:
TV: | To make sure t is in its Dir domain.
|
---|
t: | Parameter value to verify.
|
---|
Dir: | Direction. Either U or V or W.
|
---|
Returned Value:
CagdBType: TRUE if in domain, FALSE otherwise.
|
---|
Keywords:
trivariates
(triv_lib/triv_aux.c:118)
Prototype:
CagdBType TrivParamsInDomain(TrivTVStruct *TV,
CagdRType u,
CagdRType v,
CagdRType w)
Description:
Given a tri-variate and a domain - validate it.
Parameters:
TV: | To make sure (u, v, w) is in its domain.
|
---|
u, v, w: | To verify if it is in TV's parametric domain.
|
---|
Returned Value:
CagdBType: TRUE if in domain, FALSE otherwise.
|
---|
Keywords:
trivariates
(triv_lib/geomat4d.c:44)
Prototype:
int TrivPlaneFrom4Points(TrivPType Pt1,
TrivPType Pt2,
TrivPType Pt3,
TrivPType Pt4,
TrivPlaneType Plane)
Description:
Computes a hyperplane in four space through the given four points.
Based on a direct solution in Maple of:
with(linalg);
readlib(C);
d := det( matrix( [ [x - x1, y - y1, z - z1, w - w1],
[x2 - x1, y2 - y1, z2 - z1, w2 - w1],
[x3 - x2, y3 - y2, z3 - z2, w3 - w2],
[x4 - x3, y4 - y3, z4 - z3, w4 - w3]] ) );
coeff( d, x );
coeff( d, y );
coeff( d, z );
coeff( d, w );
Parameters:
Pt1, Pt2, Pt3, Pt4: | The four points the plane should go through.
|
---|
Plane: | Where the result should be placed.
|
---|
Returned Value:
int: TRUE if successful, FALSE otherwise.
|
---|
Keywords:
hyper plane
plane
(triv_lib/triveval.c:381)
Prototype:
CagdSrfStruct *TrivSrfFromMesh(TrivTVStruct *TV,
int Index,
TrivTVDirType Dir)
Description:
Extract a bivariate surface out of the given trivariate's mesh.
The provided (zero based) Index specifies which Index to extract.
Parameters:
TV: | Trivariate to extract a bivariate surface out of its mesh.
|
---|
Index: | Index of row/column/level of TV's mesh in direction Dir.
|
---|
Dir: | Direction of isosurface extraction. Either U or V or W.
|
---|
Returned Value:
CagdSrfStruct *: A bivariate surface which was extracted from TV's
Mesh. This surface is not necessarily on TV.
|
---|
Keywords:
trivariates
(triv_lib/triveval.c:202)
Prototype:
CagdSrfStruct *TrivSrfFromTV(TrivTVStruct *TV,
CagdRType t,
TrivTVDirType Dir)
Description:
Extract an isoparametric surface out of the given tensor product
trivariate.
Operations should favor the CONST_W_DIR, in which the extraction is
somewhat faster, if it is possible.
Parameters:
TV: | To extract an isoparametric surface from at parameter value t
in direction Dir.
|
---|
t: | Parameter value at which to extract the isosurface.
|
---|
Dir: | Direction of isosurface extraction. Either U or V or W.
|
---|
Returned Value:
CagdSrfStruct *: A bivariate surface which is an isosurface of TV.
|
---|
Keywords:
trivariates
(triv_lib/triveval.c:510)
Prototype:
void TrivSrfToMesh(CagdSrfStruct *Srf,
int Index,
TrivTVDirType Dir,
TrivTVStruct *TV)
Description:
Substitute a bivariate surface into a given trivariate's mesh.
The provided (zero based) Index specifies which Index to extract.
Parameters:
Srf: | Surface to substitute into the trivariate TV.
|
---|
Index: | Index of row/column/level of TV's mesh in direction Dir.
|
---|
Dir: | Direction of isosurface extraction. Either U or V or W.
|
---|
TV: | Trivariate to substitute a bivariate surface into its mesh.
|
---|
Returned Value:
Keywords:
trivariates
(triv_lib/trivmesh.c:24)
Prototype:
CagdPolylineStruct *TrivTV2CtrlMesh(TrivTVStruct *Trivar)
Description:
Extracts the control mesh of a surface as a list of polylines.
Parameters:
Srf: | To extract a control mesh from.
|
---|
Returned Value:
CagdPolylineStruct *: The control mesh of Srf.
|
---|
Keywords:
control mesh
(triv_lib/triv_aux.c:219)
Prototype:
void TrivTVBBox(TrivTVStruct *Trivar, CagdBBoxStruct *BBox)
Description:
Computes a bounding box for a trivariate freeform function.
Parameters:
Trivar: | To compute a bounding box for.
|
---|
BBox: | Where bounding information is to be saved.
|
---|
Returned Value:
Keywords:
bbox
bounding box
(triv_lib/triv_gen.c:164)
Prototype:
TrivTVStruct *TrivTVCopy(TrivTVStruct *TV)
Description:
Allocates and duplicates all slots of a trivariate structure.
Parameters:
TV: | Trivariate to duplicate
|
---|
Returned Value:
TrivTVStruct *: Duplicated trivariate.
|
---|
Keywords:
trivariates
(triv_lib/triv_gen.c:225)
Prototype:
TrivTVStruct *TrivTVCopyList(TrivTVStruct *TVList)
Description:
Duplicates a list of trivariate structures.
Parameters:
TVList: | List of trivariates to duplicate.
|
---|
Returned Value:
TrivTVStruct *: Duplicated list of trivariates.
|
---|
Keywords:
trivariates
(triv_lib/trivrais.c:35)
Prototype:
TrivTVStruct *TrivTVDegreeRaise(TrivTVStruct *TV, TrivTVDirType Dir)
Description:
Returns a new trivariate representing the same curve as TV but with its
degree raised by one.
Parameters:
TV: | To raise its degree.
|
---|
Dir: | Direction of degree raising. Either U, V or W.
|
---|
Returned Value:
TrivTVStruct *: A surface with same geometry as TV but with one
degree higher.
|
---|
Keywords:
degree raising
(triv_lib/triv_der.c:35)
Prototype:
TrivTVStruct *TrivTVDerive(TrivTVStruct *TV, TrivTVDirType Dir)
Description:
Given a trivariate, computes its partial derivative trivariate in
direction Dir.
Parameters:
TV: | Trivariate to differentiate.
|
---|
Dir: | Direction of differentiation. Either U or V or W.
|
---|
Returned Value:
TrivTVStruct *: Differentiated trivariate in direction Dir.
|
---|
Keywords:
trivariates
(triv_lib/triv_aux.c:33)
Prototype:
void TrivTVDomain(TrivTVStruct *TV,
CagdRType *UMin,
CagdRType *UMax,
CagdRType *VMin,
CagdRType *VMax,
CagdRType *WMin,
CagdRType *WMax)
Description:
Given a tri-variate, returns its parametric domain.
Parameters:
TV: | Trivariate function to consider.
|
---|
UMin, UMax: | U Domain of TV will be placed herein.
|
---|
VMin, VMax: | V Domain of TV will be placed herein.
|
---|
WMin, WMax: | W Domain of TV will be placed herein.
|
---|
Returned Value:
Keywords:
trivariates
(triv_lib/triveval.c:44)
Prototype:
CagdRType *TrivTVEval(TrivTVStruct *TV, CagdRType u, CagdRType v, CagdRType w)
Description:
Evaluates the given tensor product trivariate at a given point, by
extracting an isoparamteric surface along w and evaluating (u,v) in it.
+-----------------------+
W / /|
/ / / |
/ / U --> / | Orientation of
+-----------------------+ | control tri-variate mesh.
V | |P0 Pi-1| +
v |Pi P2i-1| /
| | /
|Pn-i Pn-1|/
+-----------------------+
Parameters:
TV: | To evaluate at given (u, v, w) parametric location.
|
---|
u, v, w: | Parametric location to evaluate TV at.
|
---|
Returned Value:
CagdRType *: A vector holding all the coefficients of all components
of the trivariate's point type. If for example trivariate
point type is P2, the W, X, and Y will be saved in the
first three locations of the returned vector. The first
location (index 0) of the returned vector is reserved for
the rational coefficient W and XYZ always starts at second
location of the returned vector (index 1).
|
---|
Keywords:
evaluation
trivariates
(triv_lib/triveval.c:163)
Prototype:
CagdRType *TrivTVEval2(TrivTVStruct *TV, CagdRType u, CagdRType v, CagdRType w)
Description:
This function is the same as TrivTVEval2 above. Cleaner, but much less
efficient.
Evaluates the given tensor product trivariate at a given point, by
extracting an isoparamteric surface along w and evaluating (u, v) in it.
Parameters:
TV: | To evaluate at given (u, v, w) parametric location.
|
---|
u, v, w: | Parametric location to evaluate TV at.
|
---|
Returned Value:
CagdRType *: A vector holding all the coefficients of all components
of the trivariate's point type. If for example trivariate
point type is P2, the W, X, and Y will be saved in the
first three locations of the returned vector. The first
location (index 0) of the returned vector is reserved for
the rational coefficient W and XYZ always starts at second
location of the returned vector (index 1).
|
---|
Keywords:
evaluation
trivariates
(triv_lib/triv_gen.c:254)
Prototype:
void TrivTVFree(TrivTVStruct *TV)
Description:
Deallocates and frees all slots of a trivariate structure.
Parameters:
Returned Value:
Keywords:
trivariates
(triv_lib/triv_gen.c:290)
Prototype:
void TrivTVFreeList(TrivTVStruct *TVList)
Description:
Deallocates and frees a list of trivariate structures.
Parameters:
TVList: | Trivariate list to free.
|
---|
Returned Value:
Keywords:
trivariates
(triv_lib/trivstrv.c:33)
Prototype:
TrivTVStruct *TrivTVFromSrfs(CagdSrfStruct *SrfList, int OtherOrder)
Description:
Constructs a trivariate using a set of surfaces. Surfaces are made to be
compatible and then each is substituted into the new trivariate's mesh as
a row.
If the OtherOrder is less than the number of curves, number of curves is
used.
A knot vector is formed with uniform open end for the other direction,
so it interpolates the first and last surfaces.
Note, however, that only the first and the last surfaces are
interpolated if OtherOrder is greater than 2.
Parameters:
SrfList: | List of surfaces to consturct a trivariate with.
|
---|
OtherOrder: | Other, third, order of trivariate.
|
---|
Returned Value:
TrivTVStruct *: Constructed trivariate from surfaces.
|
---|
Keywords:
trivar constructors
(triv_lib/triv_aux.c:246)
Prototype:
void TrivTVListBBox(TrivTVStruct *TVs, CagdBBoxStruct *BBox)
Description:
Computes a bounding box for a list of trivariate freeform function.
Parameters:
Trivars: | To compute a bounding box for.
|
---|
BBox: | Where bounding information is to be saved.
|
---|
Returned Value:
Keywords:
bbox
bounding box
(triv_lib/triv_gen.c:475)
Prototype:
void TrivTVMatTransform(TrivTVStruct *TV, CagdMType Mat)
Description:
Transforms, in place, the given TV as specified by homogeneous matrix Mat.
Parameters:
TV: | Trivariate to transform.
|
---|
Mat: | Homogeneous transformation to apply to TV.
|
---|
Returned Value:
Keywords:
trivariates
(triv_lib/triv_gen.c:35)
Prototype:
TrivTVStruct *TrivTVNew(TrivGeomType GType,
CagdPointType PType,
int ULength,
int VLength,
int WLength)
Description:
Allocates the memory required for a new trivariate.
Parameters:
GType: | Type of geometry the curve should be - Bspline, Bezier etc.
|
---|
PType: | Type of control points (E2, P3, etc.).
|
---|
ULength: | Number of control points in the U direction.
|
---|
VLength: | Number of control points in the V direction.
|
---|
WLength: | Number of control points in the W direction.
|
---|
Returned Value:
TrivTVStruct *: An uninitialized freeform trivariate.
|
---|
Keywords:
trivariates
allocation
(triv_lib/triv_ref.c:41)
Prototype:
TrivTVStruct *TrivTVRefineAtParams(TrivTVStruct *TV,
TrivTVDirType Dir,
CagdBType Replace,
CagdRType *t,
int n)
Description:
Given a trivariate, refines it at the given n knots as defined by the
vector t.
If Replace is TRUE, the values replace the current knot vector.
Returns pointer to refined TV (Note a Bezier trivariate will be
converted into a Bspline trivariate).
Parameters:
TV: | Trivariate to refine according to t in direction Dir.
|
---|
Dir: | Direction of refinement. Either U or V or W.
|
---|
Replace: | If TRUE t is a knot vector exaclt in the length of the knot
vector in direction Dir in TV and t simply replaces than knot
vector. If FALSE, the knot vector in direction Dir in TV is
refined by adding all the knots in t.
|
---|
t: | Knot vector to refine/replace the knot vector of TV in
direction Dir.
|
---|
n: | Length of vector t.
|
---|
Returned Value:
TrivTVStruct *: The refined trivariate. Always a Bspline trivariate.
|
---|
Keywords:
trivariates
(triv_lib/triv_aux.c:147)
Prototype:
TrivTVStruct *TrivTVRegionFromTV(TrivTVStruct *TV,
CagdRType t1,
CagdRType t2,
TrivTVDirType Dir)
Description:
Given a tri-variate, returns a sub-region of it.
Parameters:
TV: | To extract asub-region from.
|
---|
t1, t1: | Domain to extract from TV, in parametric direction Dir.
|
---|
Dir: | Direction to extract the sub-region. Either U or V or W.
|
---|
Returned Value:
TrivTVStruct *: A sub-region of TV from t1 to t2 in direction Dir.
|
---|
Keywords:
trivariates
(triv_lib/triv_sub.c:29)
Prototype:
TrivTVStruct *TrivTVSubdivAtParam(TrivTVStruct *TV,
CagdRType t,
TrivTVDirType Dir)
Description:
Given a tri-variate, subdivides it at parameter value t in direction
Dir.
Parameters:
TV: | Trivariate to subdivide.
|
---|
t: | Parameter to subdivide at.
|
---|
Dir: | Direction of subdivision.
|
---|
Returned Value:
TrivTVStruct *: A list of two trivariates, result of the subdivision.
|
---|
Keywords:
trivariates
(triv_lib/triv_gen.c:443)
Prototype:
void TrivTVTransform(TrivTVStruct *TV, CagdRType *Translate, CagdRType Scale)
Description:
Linearly transforms, in place, given TV as specified by Translate and
Scale.
Parameters:
TV: | Trivariate to transform.
|
---|
Translate: | Translation factor.
|
---|
Scale: | Scaling factor.
|
---|
Returned Value:
Keywords:
trivariates
(triv_lib/triv_gen.c:338)
Prototype:
TrivTriangleStruct *TrivTriangleCopy(TrivTriangleStruct *Triangle)
Description:
Allocates and duplicates all slots of a triangle structure.
Parameters:
TrivTriangleStruct *: | Triangle to duplicate.
|
---|
Returned Value:
TrivTriangleStruct *: Duplicated triangle.
|
---|
Keywords:
(triv_lib/triv_gen.c:365)
Prototype:
TrivTriangleStruct *TrivTriangleCopyList(TrivTriangleStruct *TriangleList)
Description:
Duplicates a list of triangle structures.
Parameters:
TriangleList: | List of triangle to duplicate.
|
---|
Returned Value:
TrivTriangleStruct *: Duplicated list of triangle.
|
---|
Keywords:
(triv_lib/triv_gen.c:394)
Prototype:
void TrivTriangleFree(TrivTriangleStruct *Triangle)
Description:
Deallocates and frees all slots of a triangle structure.
Parameters:
Triangle: | Triangle to free.
|
---|
Returned Value:
Keywords:
(triv_lib/triv_gen.c:416)
Prototype:
void TrivTriangleFreeList(TrivTriangleStruct *TriangleList)
Description:
Deallocates and frees a list of triangle structures.
Parameters:
TriangleList: | Triangle list to free.
|
---|
Returned Value:
Keywords:
(triv_lib/triv_gen.c:314)
Prototype:
TrivTriangleStruct *TrivTriangleNew(void)
Description:
Allocates the memory required for a new triangle.
Parameters:
None
Returned Value:
TrivTriangleStruct *: An uninitialized triangle.
|
---|
Keywords:
allocation
(triv_lib/geomat4d.c:180)
Prototype:
void TrivVectCross3Vecs(TrivVType A,
TrivVType B,
TrivVType C,
TrivVType Res)
Description:
Computes a vector in R^4 that is perpendicular to the given three vectors.
with(linalg);
readlib(C);
d := det( matrix( [ [I, J, K, L],
[A[0], A[1], A[2], A[3]],
[B[0], B[1], B[2], B[3]],
[C[0], C[1], C[2], C[3]] ] ) );
coeff( d, I );
coeff( d, J );
coeff( d, K );
coeff( d, L );
Parameters:
A, B, C: | The three vectors to compute their cross product.
|
---|
Res: | Where the output goes into.
|
---|
Returned Value:
Keywords:
cross product